1   /* Copyright 2002-2016 CS Systèmes d'Information
2    * Licensed to CS Systèmes d'Information (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.forces.drag;
18  
19  import org.hipparchus.geometry.euclidean.threed.Vector3D;
20  import org.hipparchus.util.FastMath;
21  import org.orekit.bodies.BodyShape;
22  import org.orekit.bodies.GeodeticPoint;
23  import org.orekit.errors.OrekitException;
24  import org.orekit.errors.OrekitMessages;
25  import org.orekit.frames.Frame;
26  import org.orekit.frames.Transform;
27  import org.orekit.time.AbsoluteDate;
28  import org.orekit.utils.Constants;
29  import org.orekit.utils.PVCoordinates;
30  import org.orekit.utils.PVCoordinatesProvider;
31  
32  /** This is the realization of the Jacchia-Bowman 2006 atmospheric model.
33   * <p>
34   * It is described in the paper: <br>
35   *
36   * <a href="http://sol.spacenvironment.net/~JB2006/pubs/JB2006_AIAA-6166_model.pdf">A
37   * New Empirical Thermospheric Density Model JB2006 Using New Solar Indices</a><br>
38   *
39   * <i>Bruce R. Bowman, W. Kent Tobiska and Frank A. Marcos</i> <br>
40   *
41   * AIAA 2006-6166<br>
42   *</p>
43   * <p>
44   * Two computation methods are proposed to the user:
45   * <ul>
46   * <li> one OREKIT independent and compliant with initial FORTRAN routine entry values:
47   *        {@link #getDensity(double, double, double, double, double, double, double, double, double, double, double, double, double)}. </li>
48   * <li> one compliant with OREKIT Atmosphere interface, necessary to the
49   *        {@link org.orekit.forces.drag.DragForce
50   *        drag force model} computation.</li>
51   * </ul>
52   *
53   * <p>
54   * This model provides dense output for all altitudes and positions. Output data are :
55   * <ul>
56   * <li>Exospheric Temperature above Input Position (deg K)</li>
57   * <li>Temperature at Input Position (deg K)</li>
58   * <li>Total Mass-Density at Input Position (kg/m³)</li>
59   * </ul>
60   *
61   * <p>
62   * The model needs geographical and time information to compute general values,
63   * but also needs space weather data : mean and daily solar flux, retrieved threw
64   * different indices, and planetary geomagnetic indices. <br>
65   * More information on these indices can be found on the  <a
66   * href="http://sol.spacenvironment.net/~JB2006/JB2006_index.html">
67   * official JB2006 website.</a>
68   *</p>
69   *
70   * @author Bruce R Bowman (HQ AFSPC, Space Analysis Division), Feb 2006: FORTRAN routine
71   * @author Fabien Maussion (java translation)
72   */
73  public class JB2006 implements Atmosphere {
74  
75      /** Serializable UID. */
76      private static final long serialVersionUID = -4201270765122160831L;
77  
78      /** The alpha are the thermal diffusion coefficients in equation (6). */
79      private static final double[] ALPHA = {
80          0, 0, 0, 0, 0, -0.38
81      };
82  
83      /** Natural logarithm of 10.0. */
84      private static final double AL10  = 2.3025851;
85  
86      /** Molecular weights in order: N2, O2, O, Ar, He and H. */
87      private static final double[] AMW = {
88          0,
89          28.0134, 31.9988, 15.9994, 39.9480, 4.0026, 1.00797
90      };
91  
92      /** Avogadro's number in mks units (molecules/kmol). */
93      private static final double AVOGAD = 6.02257e26;
94  
95      /** Approximate value for 2 π. */
96      private static final double TWOPI   = 6.2831853;
97  
98      /** Approximate value for π. */
99      private static final double PI      = 3.1415927;
100 
101     /** Approximate value for π / 2. */
102     private static final double PIOV2   = 1.5707963;
103 
104     /** The FRAC are the assumed sea-level volume fractions in order: N2, O2, Ar, and He. */
105     private static final double[] FRAC = {
106         0,
107         0.78110, 0.20955, 9.3400e-3, 1.2890e-5
108     };
109 
110     /** Universal gas-constant in mks units (joules/K/kmol). */
111     private static final double RSTAR = 8314.32;
112 
113     /** Value used to establish height step sizes in the regime 90km to 105km. */
114     private static final double R1 = 0.010;
115 
116     /** Value used to establish height step sizes in the regime 105km to 500km. */
117     private static final double R2 = 0.025;
118 
119     /** Value used to establish height step sizes in the regime above 500km. */
120     private static final double R3 = 0.075;
121 
122     /** Weights for the Newton-Cotes five-points quadrature formula. */
123     private static final double[] WT = {
124         0,
125         0.311111111111111, 1.422222222222222,
126         0.533333333333333, 1.422222222222222,
127         0.311111111111111
128     };
129 
130     /** Coefficients for high altitude density correction. */
131     private static final double[] CHT = {
132         0, 0.22, -0.20e-02, 0.115e-02, -0.211e-05
133     };
134 
135     /** FZ global model values (1978-2004 fit).  */
136     private static final double[] FZM = {
137         0,
138         0.111613e+00, -0.159000e-02, 0.126190e-01,
139         -0.100064e-01, -0.237509e-04, 0.260759e-04
140     };
141 
142     /** GT global model values (1978-2004 fit). */
143     private static final double[] GTM = {
144         0,
145         -0.833646e+00, -0.265450e+00, 0.467603e+00, -0.299906e+00,
146         -0.105451e+00, -0.165537e-01, -0.380037e-01, -0.150991e-01,
147         -0.541280e-01, 0.119554e-01, 0.437544e-02, -0.369016e-02,
148         0.206763e-02, -0.142888e-02, -0.867124e-05, 0.189032e-04,
149         0.156988e-03,  0.491286e-03, -0.391484e-04, -0.126854e-04,
150         0.134078e-04, -0.614176e-05, 0.343423e-05
151     };
152 
153     /** XAMBAR relative data. */
154     private static final double[] CXAMB = {
155         0,
156         28.15204, -8.5586e-2, +1.2840e-4, -1.0056e-5,
157         -1.0210e-5, +1.5044e-6, +9.9826e-8
158     };
159 
160     /** DTSUB relative data. */
161     private static final double[] BDT_SUB = {
162         0,
163         -0.457512297e+01, -0.512114909e+01, -0.693003609e+02,
164         0.203716701e+03,  0.703316291e+03, -0.194349234e+04,
165         0.110651308e+04, -0.174378996e+03,  0.188594601e+04,
166         -0.709371517e+04,  0.922454523e+04, -0.384508073e+04,
167         -0.645841789e+01,  0.409703319e+02, -0.482006560e+03,
168         0.181870931e+04, -0.237389204e+04,  0.996703815e+03,
169         0.361416936e+02
170     };
171 
172     /** DTSUB relative data. */
173     private static final double[] CDT_SUB = {
174         0,
175         -0.155986211e+02, -0.512114909e+01, -0.693003609e+02,
176         0.203716701e+03,  0.703316291e+03, -0.194349234e+04,
177         0.110651308e+04, -0.220835117e+03,  0.143256989e+04,
178         -0.318481844e+04,  0.328981513e+04, -0.135332119e+04,
179         0.199956489e+02, -0.127093998e+02,  0.212825156e+02,
180         -0.275555432e+01,  0.110234982e+02,  0.148881951e+03,
181         -0.751640284e+03,  0.637876542e+03,  0.127093998e+02,
182         -0.212825156e+02,  0.275555432e+01
183     };
184 
185     /** Temperatures.
186      *  <ul>
187      *  <li>TEMP(1): Exospheric Temperature above Input Position (deg K)</li>
188      *  <li>TEMP(2): Temperature at Input Position (deg K)</li>
189      *  </ul>
190      */
191     private double[] temp = new double[3];
192 
193     /** Total Mass-Density at Input Position (kg/m³). */
194     private double rho;
195 
196     /** Sun position. */
197     private PVCoordinatesProvider sun;
198 
199     /** External data container. */
200     private JB2006InputParameters inputParams;
201 
202     /** Earth body shape. */
203     private BodyShape earth;
204 
205     /** Constructor with space environment information for internal computation.
206      * @param parameters the solar and magnetic activity data
207      * @param sun the sun position
208      * @param earth the earth body shape
209      */
210     public JB2006(final JB2006InputParameters parameters,
211                   final PVCoordinatesProvider sun, final BodyShape earth) {
212         this.earth = earth;
213         this.sun = sun;
214         this.inputParams = parameters;
215     }
216 
217     /** {@inheritDoc} */
218     public Frame getFrame() {
219         return earth.getBodyFrame();
220     }
221 
222     /** Get the local density with initial entries.
223      * @param dateMJD date and time, in modified julian days and fraction
224      * @param sunRA Right Ascension of Sun (radians)
225      * @param sunDecli Declination of Sun (radians)
226      * @param satLon Right Ascension of position (radians)
227      * @param satLat Geocentric latitude of position (radians)
228      * @param satAlt Height of position (m)
229      * @param f10 10.7-cm Solar flux (1e<sup>-22</sup>*Watt/(m²*Hertz)).
230      *            Tabular time 1.0 day earlier
231      * @param f10B 10.7-cm Solar Flux, averaged 81-day centered on the input time
232      * @param ap Geomagnetic planetary 3-hour index A<sub>p</sub>
233      *            for a tabular time 6.7 hours earlier
234      * @param s10 EUV index (26-34 nm) scaled to F10. Tabular time 1 day earlier.
235      * @param s10B UV 81-day averaged centered index
236      * @param xm10 MG2 index scaled to F10
237      * @param xm10B MG2 81-day ave. centered index. Tabular time 5.0 days earlier.
238      * @return total mass-Density at input position (kg/m³)
239      * @exception OrekitException if altitude is below 90 km
240      */
241     public double getDensity(final double dateMJD, final double sunRA, final double sunDecli,
242                              final double satLon, final double satLat, final double satAlt,
243                              final double f10, final double f10B, final double ap,
244                              final double s10, final double s10B, final double xm10, final double xm10B)
245         throws OrekitException {
246 
247         if (satAlt < 90000) {
248             throw new OrekitException(OrekitMessages.ALTITUDE_BELOW_ALLOWED_THRESHOLD,
249                                       satAlt, 90000.0);
250         }
251         final double scaledSatAlt = satAlt / 1000.0;
252 
253         // Equation (14)
254         final double tc = 379 + 3.353 * f10B + 0.358 * (f10 - f10B) +
255                           2.094 * (s10 - s10B) + 0.343 * (xm10 - xm10B);
256 
257         // Equation (15)
258         final double eta =   0.5 * FastMath.abs(satLat - sunDecli);
259         final double theta = 0.5 * FastMath.abs(satLat + sunDecli);
260 
261         // Equation (16)
262         final double h     = satLon - sunRA;
263         final double tau   = h - 0.64577182 + 0.10471976 * FastMath.sin(h + 0.75049158);
264         double solTimeHour = FastMath.toDegrees(h + PI) / 15.0;
265         if (solTimeHour >= 24) {
266             solTimeHour = solTimeHour - 24.;
267         }
268         if (solTimeHour < 0) {
269             solTimeHour = solTimeHour + 24.;
270         }
271 
272         // Equation (17)
273         final double C = FastMath.pow(FastMath.cos(eta), 2.5);
274         final double S = FastMath.pow(FastMath.sin(theta), 2.5);
275         final double tmp = FastMath.abs(FastMath.cos(0.5 * tau));
276         final double DF = S + (C - S) * tmp * tmp * tmp;
277         final double TSUBL = tc * (1. + 0.31 * DF);
278 
279         // Equation (18)
280         final double EXPAP = FastMath.exp(-0.08 * ap);
281         final double DTG = ap + 100. * (1. - EXPAP);
282 
283         // Compute correction to dTc for local solar time and lat correction
284         final double DTCLST = dTc(f10, solTimeHour, satLat, scaledSatAlt);
285 
286         // Compute the local exospheric temperature.
287         final double TINF = TSUBL + DTG + DTCLST;
288         temp[1] = TINF;
289 
290         // Equation (9)
291         final double TSUBX = 444.3807 + 0.02385 * TINF - 392.8292 * FastMath.exp(-0.0021357 * TINF);
292 
293         // Equation (11)
294         final double GSUBX = 0.054285714 * (TSUBX - 183.);
295 
296         // The TC array will be an argument in the call to
297         // XLOCAL, which evaluates Equation (10) or Equation (13)
298         final double[] TC = new double[4];
299         TC[0] = TSUBX;
300         TC[1] = GSUBX;
301 
302         //   A AND GSUBX/A OF Equation (13)
303         TC[2] = (TINF - TSUBX) / PIOV2;
304         TC[3] = GSUBX / TC[2];
305 
306         // Equation (5)
307         final double Z1 = 90.;
308         final double Z2 = FastMath.min(scaledSatAlt, 105.0);
309         double AL = FastMath.log(Z2 / Z1);
310         int N = (int) FastMath.floor(AL / R1) + 1;
311         double ZR = FastMath.exp(AL / N);
312         final double AMBAR1 = xAmbar(Z1);
313         final double TLOC1 = xLocal(Z1, TC);
314         double ZEND   = Z1;
315         double SUM2   = 0.;
316         double AIN    = AMBAR1 * xGrav(Z1) / TLOC1;
317         double AMBAR2 = 0;
318         double TLOC2  = 0;
319         double Z      = 0;
320         double GRAVL  = 0;
321 
322         for (int i = 1; i <= N; ++i) {
323             Z = ZEND;
324             ZEND = ZR * Z;
325             final double DZ = 0.25 * (ZEND - Z);
326             double SUM1 = WT[1] * AIN;
327             for (int j = 2; j <= 5; ++j) {
328                 Z += DZ;
329                 AMBAR2 = xAmbar(Z);
330                 TLOC2  = xLocal(Z, TC);
331                 GRAVL  = xGrav(Z);
332                 AIN    = AMBAR2 * GRAVL / TLOC2;
333                 SUM1  += WT[j] * AIN;
334             }
335             SUM2 = SUM2 + DZ * SUM1;
336         }
337         final double FACT1 = 1000.0 / RSTAR;
338         rho = 3.46e-6 * AMBAR2 * TLOC1 * FastMath.exp(-FACT1 * SUM2) / (AMBAR1 * TLOC2);
339 
340         // Equation (2)
341         final double ANM = AVOGAD * rho;
342         double AN  = ANM / AMBAR2;
343 
344         // Equation (3)
345         double FACT2  = ANM / 28.960;
346         final double[] ALN = new double[7];
347         ALN[1] = FastMath.log(FRAC[1] * FACT2);
348         ALN[4] = FastMath.log(FRAC[3] * FACT2);
349         ALN[5] = FastMath.log(FRAC[4] * FACT2);
350 
351         // Equation (4)
352         ALN[2] = FastMath.log(FACT2 * (1. + FRAC[2]) - AN);
353         ALN[3] = FastMath.log(2. * (AN - FACT2));
354 
355         if (scaledSatAlt <= 105.0) {
356             temp[2] = TLOC2;
357             // Put in negligible hydrogen for use in DO-LOOP 13
358             ALN[6] = ALN[5] - 25.0;
359         } else {
360             // Equation (6)
361             final double Z3 = FastMath.min(scaledSatAlt, 500.0);
362             AL   = FastMath.log(Z3 / Z);
363             N    = (int) FastMath.floor(AL / R2) + 1;
364             ZR   = FastMath.exp(AL / N);
365             SUM2 = 0.;
366             AIN  = GRAVL / TLOC2;
367 
368             double TLOC3 = 0;
369             for (int I = 1; I <= N; ++I) {
370                 Z = ZEND;
371                 ZEND = ZR * Z;
372                 final double DZ = 0.25 * (ZEND - Z);
373                 double SUM1 = WT[1] * AIN;
374                 for (int J = 2; J <= 5; ++J) {
375                     Z    += DZ;
376                     TLOC3 = xLocal(Z, TC);
377                     GRAVL = xGrav(Z);
378                     AIN   = GRAVL / TLOC3;
379                     SUM1  = SUM1 + WT[J] * AIN;
380                 }
381                 SUM2 = SUM2 + DZ * SUM1;
382             }
383 
384             final double Z4 = FastMath.max(scaledSatAlt, 500.0);
385             AL = FastMath.log(Z4 / Z);
386             double R = R2;
387             if (scaledSatAlt > 500.0) {
388                 R = R3;
389             }
390             N = (int) FastMath.floor(AL / R) + 1;
391             ZR = FastMath.exp(AL / N);
392             double SUM3 = 0.;
393             double TLOC4 = 0;
394             for (int I = 1; I <= N; ++I) {
395                 Z = ZEND;
396                 ZEND = ZR * Z;
397                 final double DZ = 0.25 * (ZEND - Z);
398                 double SUM1 = WT[1] * AIN;
399                 for (int J = 2; J <= 5; ++J) {
400                     Z    += DZ;
401                     TLOC4 = xLocal(Z, TC);
402                     GRAVL = xGrav(Z);
403                     AIN   = GRAVL / TLOC4;
404                     SUM1  = SUM1 + WT[J] * AIN;
405                 }
406                 SUM3 = SUM3 + DZ * SUM1;
407             }
408             final double ALTR;
409             final double HSIGN;
410             if (scaledSatAlt <= 500.) {
411                 temp[2] = TLOC3;
412                 ALTR = FastMath.log(TLOC3 / TLOC2);
413                 FACT2 = FACT1 * SUM2;
414                 HSIGN = 1.0;
415 
416             } else {
417                 temp[2] = TLOC4;
418                 ALTR = FastMath.log(TLOC4 / TLOC2);
419                 FACT2 = FACT1 * (SUM2 + SUM3);
420                 HSIGN = -1.0;
421             }
422             for (int I = 1; I <= 5; ++I) {
423                 ALN[I] = ALN[I] - (1.0 + ALPHA[I]) * ALTR - FACT2 * AMW[I];
424             }
425 
426             // Equation (7) - Note that in CIRA72, AL10T5 = DLOG10(T500)
427             final double AL10T5 = FastMath.log(TINF) / FastMath.log(10);
428             final double ALNH5 = (5.5 * AL10T5 - 39.40) * AL10T5 + 73.13;
429             ALN[6] = AL10 * (ALNH5 + 6.) + HSIGN * (FastMath.log(TLOC4 / TLOC3) + FACT1 * SUM3 * AMW[6]);
430 
431         }
432 
433         // Equation (24)  - J70 Seasonal-Latitudinal Variation
434         final double CAPPHI = ((dateMJD - 36204.0) / 365.2422) % 1;
435         final int signum = (satLat >= 0) ? 1 : -1;
436         final double sinLat = FastMath.sin(satLat);
437         final double DLRSL = 0.02 * (scaledSatAlt - 90.) * FastMath.exp(-0.045 * (scaledSatAlt - 90.)) *
438                              signum * FastMath.sin(TWOPI * CAPPHI + 1.72) * sinLat * sinLat;
439 
440         // Equation (23) - Computes the semiannual variation
441         double DLRSA = 0;
442         if (Z < 2000.0) {
443             final double D1950 = dateMJD - 33281.0;
444             // Use new semiannual model DELTA LOG RHO
445             DLRSA = semian(dayOfYear(D1950), scaledSatAlt, f10B);
446         }
447 
448         // Sum the delta-log-rhos and apply to the number densities.
449         // In CIRA72 the following equation contains an actual sum,
450         // namely DLR = AL10 * (DLRGM + DLRSA + DLRSL)
451         // However, for Jacchia 70, there is no DLRGM or DLRSA.
452         final double DLR = AL10 * (DLRSL + DLRSA);
453         for (int i = 1; i <= 6; ++i) {
454             ALN[i] += DLR;
455         }
456 
457         // Compute mass-density and mean-molecular-weight and
458         // convert number density logs from natural to common.
459 
460         double SUMNM = 0.0;
461 
462         for (int I = 1; I <= 6; ++I) {
463             AN = FastMath.exp(ALN[I]);
464             SUMNM += AN * AMW[I];
465         }
466 
467         rho = SUMNM / AVOGAD;
468 
469         // Compute the high altitude exospheric density correction factor
470         double FEX = 1.0;
471         if ((scaledSatAlt >= 1000.0) && (scaledSatAlt < 1500.0)) {
472             final double ZETA   = (scaledSatAlt - 1000.) * 0.002;
473             final double ZETA2  =  ZETA * ZETA;
474             final double ZETA3  =  ZETA * ZETA2;
475             final double F15C   = CHT[1] + CHT[2] * f10B + CHT[3] * 1500.0 + CHT[4] * f10B * 1500.0;
476             final double F15C_ZETA = (CHT[3] + CHT[4] * f10B) * 500.0;
477             final double FEX2   = 3.0 * F15C - F15C_ZETA - 3.0;
478             final double FEX3   = F15C_ZETA - 2.0 * F15C + 2.0;
479             FEX    = 1.0 + FEX2 * ZETA2 + FEX3 * ZETA3;
480         }
481         if (scaledSatAlt >= 1500.0) {
482             FEX    = CHT[1] + CHT[2] * f10B + CHT[3] * scaledSatAlt + CHT[4] * f10B * scaledSatAlt;
483         }
484 
485         // Apply the exospheric density correction factor.
486         rho  *= FEX;
487 
488         return rho;
489 
490     }
491 
492     /** Compute daily temperature correction for Jacchia-Bowman model.
493      * @param f10 solar flux index
494      * @param solTimeHour local solar time (hours 0-23.999)
495      * @param satLat sat lat (radians)
496      * @param satAlt height (km)
497      * @return dTc correction
498      */
499     private static double dTc(final double f10, final double solTimeHour,
500                               final double satLat, final double satAlt) {
501 
502         double dTc = 0;
503         final double tx  = solTimeHour / 24.0;
504         final double tx2 = tx * tx;
505         final double tx3 = tx2 * tx;
506         final double tx4 = tx3 * tx;
507         final double tx5 = tx4 * tx;
508         final double ycs = FastMath.cos(satLat);
509         final double f   = (f10 - 100.0) / 100.0;
510         double h;
511         double sum;
512 
513         // Calculates dTc
514         if ((satAlt >= 120) && (satAlt <= 200)) {
515             final double DTC200 = CDT_SUB[17] + CDT_SUB[18] * tx * ycs + CDT_SUB[19] * tx2 * ycs +
516                                   CDT_SUB[20] * tx3 * ycs + CDT_SUB[21] * f * ycs + CDT_SUB[22] * tx * f * ycs +
517                                   CDT_SUB[23] * tx2 * f * ycs;
518             sum = CDT_SUB[1] + BDT_SUB[2] * f + CDT_SUB[3] * tx * f + CDT_SUB[4] * tx2 * f +
519                   CDT_SUB[5] * tx3 * f + CDT_SUB[6] * tx4 * f + CDT_SUB[7] * tx5 * f +
520                   CDT_SUB[8] * tx * ycs + CDT_SUB[9] * tx2 * ycs + CDT_SUB[10] * tx3 * ycs +
521                   CDT_SUB[11] * tx4 * ycs + CDT_SUB[12] * tx5 * ycs + CDT_SUB[13] * ycs +
522                   CDT_SUB[14] * f * ycs + CDT_SUB[15] * tx * f * ycs  + CDT_SUB[16] * tx2 * f * ycs;
523             final double DTC200DZ = sum;
524             final double CC  = 3.0 * DTC200 - DTC200DZ;
525             final double DD  = DTC200 - CC;
526             final double ZP  = (satAlt - 120.0) / 80.0;
527             dTc = CC * ZP * ZP + DD * ZP * ZP * ZP;
528         }
529 
530         if ((satAlt > 200.0) && (satAlt <= 240.0)) {
531             h = (satAlt - 200.0) / 50.0;
532             sum = CDT_SUB[1] * h + BDT_SUB[2] * f * h + CDT_SUB[3] * tx * f * h + CDT_SUB[4] * tx2 * f * h +
533                   CDT_SUB[5] * tx3 * f * h + CDT_SUB[6] * tx4 * f * h + CDT_SUB[7] * tx5 * f * h +
534                   CDT_SUB[8] * tx * ycs * h + CDT_SUB[9] * tx2 * ycs * h + CDT_SUB[10] * tx3 * ycs * h +
535                   CDT_SUB[11] * tx4 * ycs * h + CDT_SUB[12] * tx5 * ycs * h + CDT_SUB[13] * ycs * h +
536                   CDT_SUB[14] * f * ycs * h + CDT_SUB[15] * tx * f * ycs * h  + CDT_SUB[16] * tx2 * f * ycs * h +
537                   CDT_SUB[17] + CDT_SUB[18] * tx * ycs + CDT_SUB[19] * tx2 * ycs +
538                   CDT_SUB[20] * tx3 * ycs + CDT_SUB[21] * f * ycs + CDT_SUB[22] * tx * f * ycs +
539                   CDT_SUB[23] * tx2 * f * ycs;
540             dTc = sum;
541         }
542 
543         if ((satAlt > 240.0) && (satAlt <= 300.0)) {
544             h = 40.0 / 50.0;
545             sum = CDT_SUB[1] * h + BDT_SUB[2] * f * h + CDT_SUB[3] * tx * f * h + CDT_SUB[4] * tx2 * f * h +
546                   CDT_SUB[5] * tx3 * f * h + CDT_SUB[6] * tx4 * f * h + CDT_SUB[7] * tx5 * f * h +
547                   CDT_SUB[8] * tx * ycs * h + CDT_SUB[9] * tx2 * ycs * h + CDT_SUB[10] * tx3 * ycs * h +
548                   CDT_SUB[11] * tx4 * ycs * h + CDT_SUB[12] * tx5 * ycs * h + CDT_SUB[13] * ycs * h +
549                   CDT_SUB[14] * f * ycs * h + CDT_SUB[15] * tx * f * ycs * h  + CDT_SUB[16] * tx2 * f * ycs * h +
550                   CDT_SUB[17] + CDT_SUB[18] * tx * ycs + CDT_SUB[19] * tx2 * ycs +
551                   CDT_SUB[20] * tx3 * ycs + CDT_SUB[21] * f * ycs + CDT_SUB[22] * tx * f * ycs +
552                   CDT_SUB[23] * tx2 * f * ycs;
553             final double AA = sum;
554             final double BB = CDT_SUB[1] + BDT_SUB[2] * f + CDT_SUB[3] * tx * f + CDT_SUB[4] * tx2 * f +
555                         CDT_SUB[5] * tx3 * f + CDT_SUB[6] * tx4 * f + CDT_SUB[7] * tx5 * f +
556                         CDT_SUB[8] * tx * ycs + CDT_SUB[9] * tx2 * ycs + CDT_SUB[10] * tx3 * ycs +
557                         CDT_SUB[11] * tx4 * ycs + CDT_SUB[12] * tx5 * ycs + CDT_SUB[13] * ycs +
558                         CDT_SUB[14] * f * ycs + CDT_SUB[15] * tx * f * ycs + CDT_SUB[16] * tx2 * f * ycs;
559             h   = 300.0 / 100.0;
560             sum = BDT_SUB[1] + BDT_SUB[2] * f  + BDT_SUB[3] * tx * f         + BDT_SUB[4] * tx2 * f +
561                   BDT_SUB[5] * tx3 * f      + BDT_SUB[6] * tx4 * f      + BDT_SUB[7] * tx5 * f +
562                   BDT_SUB[8] * tx * ycs       + BDT_SUB[9] * tx2 * ycs    + BDT_SUB[10] * tx3 * ycs +
563                   BDT_SUB[11] * tx4 * ycs   + BDT_SUB[12] * tx5 * ycs   + BDT_SUB[13] * h * ycs +
564                   BDT_SUB[14] * tx * h * ycs    + BDT_SUB[15] * tx2 * h * ycs + BDT_SUB[16] * tx3 * h * ycs +
565                   BDT_SUB[17] * tx4 * h * ycs + BDT_SUB[18] * tx5 * h * ycs + BDT_SUB[19] * ycs;
566             final double DTC300 = sum;
567             sum = BDT_SUB[13] * ycs +
568                   BDT_SUB[14] * tx * ycs + BDT_SUB[15] * tx2 * ycs + BDT_SUB[16] * tx3 * ycs +
569                   BDT_SUB[17] * tx4 * ycs + BDT_SUB[18] * tx5 * ycs;
570             final double DTC300DZ = sum;
571             final double CC = 3.0 * DTC300 - DTC300DZ - 3.0 * AA - 2.0 * BB;
572             final double  DD = DTC300 - AA - BB - CC;
573             final double ZP  = (satAlt - 240.0) / 60.0;
574             dTc = AA + BB * ZP + CC * ZP * ZP + DD * ZP * ZP * ZP;
575         }
576 
577         if ((satAlt > 300.0) && (satAlt <= 600.0)) {
578             h   = satAlt / 100.0;
579             sum = BDT_SUB[1]    + BDT_SUB[2] * f  + BDT_SUB[3] * tx * f         + BDT_SUB[4] * tx2 * f +
580                   BDT_SUB[5] * tx3 * f      + BDT_SUB[6] * tx4 * f      + BDT_SUB[7] * tx5 * f +
581                   BDT_SUB[8] * tx * ycs       + BDT_SUB[9] * tx2 * ycs    + BDT_SUB[10] * tx3 * ycs +
582                   BDT_SUB[11] * tx4 * ycs   + BDT_SUB[12] * tx5 * ycs   + BDT_SUB[13] * h * ycs +
583                   BDT_SUB[14] * tx * h * ycs    + BDT_SUB[15] * tx2 * h * ycs + BDT_SUB[16] * tx3 * h * ycs +
584                   BDT_SUB[17] * tx4 * h * ycs + BDT_SUB[18] * tx5 * h * ycs + BDT_SUB[19] * ycs;
585             dTc = sum;
586         }
587 
588         if ((satAlt > 600.0) && (satAlt <= 800.0)) {
589             final double ZP = (satAlt - 600.0) / 100.0;
590             final double HP = 600.0 / 100.0;
591             final double AA  = BDT_SUB[1]    + BDT_SUB[2] * f  + BDT_SUB[3] * tx * f         + BDT_SUB[4] * tx2 * f +
592                                BDT_SUB[5] * tx3 * f      + BDT_SUB[6] * tx4 * f      + BDT_SUB[7] * tx5 * f +
593                                BDT_SUB[8] * tx * ycs       + BDT_SUB[9] * tx2 * ycs    + BDT_SUB[10] * tx3 * ycs +
594                                BDT_SUB[11] * tx4 * ycs   + BDT_SUB[12] * tx5 * ycs   + BDT_SUB[13] * HP * ycs +
595                                BDT_SUB[14] * tx * HP * ycs   + BDT_SUB[15] * tx2 * HP * ycs + BDT_SUB[16] * tx3 * HP * ycs +
596                                BDT_SUB[17] * tx4 * HP * ycs + BDT_SUB[18] * tx5 * HP * ycs + BDT_SUB[19] * ycs;
597             final double BB  = BDT_SUB[13] * ycs +
598                                BDT_SUB[14] * tx * ycs    + BDT_SUB[15] * tx2 * ycs + BDT_SUB[16] * tx3 * ycs +
599                                BDT_SUB[17] * tx4 * ycs + BDT_SUB[18] * tx5 * ycs;
600             final double CC  = -(3.0 * AA + 4.0 * BB) / 4.0;
601             final double DD  = (AA + BB) / 4.0;
602             dTc = AA + BB * ZP + CC * ZP * ZP + DD * ZP * ZP * ZP;
603         }
604 
605         return dTc;
606     }
607 
608     /** Evaluates Equation (1).
609      * @param z altitude
610      * @return equation (1) value
611      */
612     private static double xAmbar(final double z) {
613         final double dz = z - 100.;
614         double amb = CXAMB[7];
615         for (int i = 6; i >= 1; --i) {
616             amb = dz * amb + CXAMB[i];
617         }
618         return amb;
619     }
620 
621     /**  Evaluates Equation (10) or Equation (13), depending on Z.
622      * @param z altitude
623      * @param TC tc array ???
624      * @return equation (10) value
625      */
626     private static double xLocal(final double z, final double[] TC) {
627         final double dz = z - 125;
628         if (dz <= 0) {
629             return ((-9.8204695e-6 * dz - 7.3039742e-4) * dz * dz + 1.0) * dz * TC[1] + TC[0];
630         } else {
631             return TC[0] + TC[2] * FastMath.atan(TC[3] * dz * (1 + 4.5e-6 * FastMath.pow(dz, 2.5)));
632         }
633     }
634 
635     /** Evaluates Equation (8) of gravity field.
636      * @param z altitude
637      * @return the gravity field
638      */
639     private static double xGrav(final double z) {
640         final double temp = 1.0 + z / 6356.766;
641         return 9.80665 / (temp * temp);
642     }
643 
644     /** Compute semi-annual variation (delta log(rho)).
645      * @param day day of year
646      * @param height height (km)
647      * @param f10Bar average 81-day centered f10
648      * @return semi-annual variation
649      */
650     private static double semian (final double day, final double height, final double f10Bar) {
651 
652         final double f10Bar2 = f10Bar * f10Bar;
653         final double htz = height / 1000.0;
654 
655         // SEMIANNUAL AMPLITUDE
656         final double fzz = FZM[1] + FZM[2] * f10Bar  + FZM[3] * f10Bar * htz +
657                            FZM[4] * f10Bar * htz * htz + FZM[5] * f10Bar * f10Bar * htz +
658                            FZM[6] * f10Bar * f10Bar * htz * htz;
659 
660         // SEMIANNUAL PHASE FUNCTION
661         final double tau   = TWOPI * (day - 1.0) / 365;
662         final double sin1P = FastMath.sin(tau);
663         final double cos1P = FastMath.cos(tau);
664         final double sin2P = FastMath.sin(2.0 * tau);
665         final double cos2P = FastMath.cos(2.0 * tau);
666         final double sin3P = FastMath.sin(3.0 * tau);
667         final double cos3P = FastMath.cos(3.0 * tau);
668         final double sin4P = FastMath.sin(4.0 * tau);
669         final double cos4P = FastMath.cos(4.0 * tau);
670         final double gtz = GTM[1] + GTM[2] * sin1P + GTM[3] * cos1P +
671                            GTM[4] * sin2P + GTM[5] * cos2P +
672                            GTM[6] * sin3P + GTM[7] * cos3P +
673                            GTM[8] * sin4P + GTM[9] * cos4P +
674                            GTM[10] * f10Bar + GTM[11] * f10Bar * sin1P + GTM[12] * f10Bar * cos1P +
675                            GTM[13] * f10Bar * sin2P + GTM[14] * f10Bar * cos2P +
676                            GTM[15] * f10Bar * sin3P + GTM[16] * f10Bar * cos3P +
677                            GTM[17] * f10Bar * sin4P + GTM[18] * f10Bar * cos4P +
678                            GTM[19] * f10Bar2  + GTM[20] * f10Bar2  * sin1P + GTM[21] * f10Bar2  * cos1P +
679                            GTM[22] * f10Bar2  * sin2P + GTM[23] * f10Bar2  * cos2P;
680 
681         return FastMath.max(1.0e-6, fzz) * gtz;
682 
683     }
684 
685     /** Compute day of year.
686      * @param d1950 (days since 1950)
687      * @return the number days in year
688      */
689     private static double dayOfYear(final double d1950) {
690 
691         int iyday = (int) d1950;
692         final double frac = d1950 - iyday;
693         iyday = iyday + 364;
694 
695         int itemp = iyday / 1461;
696 
697         iyday = iyday - itemp * 1461;
698         itemp = iyday / 365;
699         if (itemp >= 3) {
700             itemp = 3;
701         }
702         iyday = iyday - 365 * itemp + 1;
703         return iyday + frac;
704     }
705 
706     // OUTPUT:
707 
708     /** Get the exospheric temperature above input position.
709      * {@link #getDensity(double, double, double, double, double, double, double, double, double, double, double, double, double)}
710      * <b> must </b> must be called before calling this function.
711      * @return the exospheric temperature (deg K)
712      */
713     public double getExosphericTemp() {
714         return temp[1];
715     }
716 
717     /** Get the temperature at input position.
718      * {@link #getDensity(double, double, double, double, double, double, double, double, double, double, double, double, double)}
719      * <b> must </b> must be called before calling this function.
720      * @return the local temperature (deg K)
721      */
722     public double getLocalTemp() {
723         return temp[2];
724     }
725 
726     /** Get the local density.
727      * @param date current date
728      * @param position current position in frame
729      * @param frame the frame in which is defined the position
730      * @return local density (kg/m³)
731      * @exception OrekitException if date is out of range of solar activity
732      */
733     public double getDensity(final AbsoluteDate date, final Vector3D position,
734                              final Frame frame)
735         throws OrekitException {
736         // check if data are available :
737         if (date.compareTo(inputParams.getMaxDate()) > 0 ||
738             date.compareTo(inputParams.getMinDate()) < 0) {
739             throw new OrekitException(OrekitMessages.NO_SOLAR_ACTIVITY_AT_DATE,
740                                       date, inputParams.getMinDate(), inputParams.getMaxDate());
741         }
742 
743         // compute modified julian days date
744         final double dateMJD = date.durationFrom(AbsoluteDate.MODIFIED_JULIAN_EPOCH) / Constants.JULIAN_DAY;
745 
746         // compute geodetic position
747         final GeodeticPoint inBody = earth.transform(position, frame, date);
748 
749         // compute sun position
750         final Frame ecef = earth.getBodyFrame();
751         final GeodeticPoint sunInBody =
752             earth.transform(sun.getPVCoordinates(date, ecef).getPosition(), ecef, date);
753         return getDensity(dateMJD,
754                           sunInBody.getLongitude(), sunInBody.getLatitude(),
755                           inBody.getLongitude(), inBody.getLatitude(),
756                           inBody.getAltitude(), inputParams.getF10(date),
757                           inputParams.getF10B(date),
758                           inputParams.getAp(date), inputParams.getS10(date),
759                           inputParams.getS10B(date), inputParams.getXM10(date),
760                           inputParams.getXM10B(date));
761     }
762 
763     /** Get the inertial velocity of atmosphere molecules.
764      * Here the case is simplified : atmosphere is supposed to have a null velocity
765      * in earth frame.
766      * @param date current date
767      * @param position current position in frame
768      * @param frame the frame in which is defined the position
769      * @return velocity (m/s) (defined in the same frame as the position)
770      * @exception OrekitException if some frame conversion cannot be performed
771      */
772     public Vector3D getVelocity(final AbsoluteDate date, final Vector3D position,
773                                 final Frame frame)
774         throws OrekitException {
775         final Transform bodyToFrame = earth.getBodyFrame().getTransformTo(frame, date);
776         final Vector3D posInBody = bodyToFrame.getInverse().transformPosition(position);
777         final PVCoordinates pvBody = new PVCoordinates(posInBody, new Vector3D(0, 0, 0));
778         final PVCoordinates pvFrame = bodyToFrame.transformPVCoordinates(pvBody);
779         return pvFrame.getVelocity();
780     }
781 
782 }